predict_2lsom Subroutine

public subroutine predict_2lsom(kohonen_map, input_data, map_output)

Subroutine to make a prediction from a trained two_level self_organized_map

Type Bound

two_level_self_organizing_map

Arguments

Type IntentOptional Attributes Name
class(two_level_self_organizing_map) :: kohonen_map

A two_level_self_organizing_map object

type(kohonen_pattern), intent(inout), dimension(:) :: input_data

A kohonen_pattern array

integer, intent(out), dimension(:,:) :: map_output

An integer array


Calls

proc~~predict_2lsom~~CallsGraph proc~predict_2lsom two_level_self_organizing_map%predict_2lsom float float proc~predict_2lsom->float none~distance~8 kohonen_prototype%distance proc~predict_2lsom->none~distance~8 none~get_prototype kohonen_prototype%get_prototype proc~predict_2lsom->none~get_prototype proc~kohonen_pattern_accessor kohonen_pattern%kohonen_pattern_accessor proc~predict_2lsom->proc~kohonen_pattern_accessor none~distance~8->none~get_prototype calculate calculate none~distance~8->calculate

Variables

Type Visibility Attributes Name Initial
integer, public :: ipattern
integer, public :: ihit
integer, public :: jhit
integer, public :: khit
integer, public :: ix
integer, public :: iy
integer, public :: iz
integer, public :: number_variables
integer, public :: ic
real(kind=wp), public :: dist_hit
real(kind=wp), public :: dist
type(kohonen_prototype), public :: current_prototype
type(kohonen_prototype), public :: predict_grid_prototype
type(kohonen_prototype), public :: predict_cluster_prototype
real(kind=wp), public, dimension(kohonen_map%number_variables1,kohonen_map%number_variables2) :: current_values
real(kind=wp), public, dimension(kohonen_map%number_clusters) :: distance_units

Source Code

   subroutine predict_2lsom(kohonen_map,input_data,map_output)
   !========================================================================================
!! Subroutine to make a prediction from a trained two_level self_organized_map 
      class(two_level_self_organizing_map) :: kohonen_map
!! A `two_level_self_organizing_map` object
      type(kohonen_pattern),dimension(:),intent(inout) :: input_data
!! A `kohonen_pattern` array
      integer,dimension(:,:),intent(out) :: map_output
!! An integer array
      integer :: ipattern,ihit,jhit,khit,ix,iy,iz,number_variables,ic
      real(kind=wp) :: dist_hit,dist
      type(kohonen_prototype) :: current_prototype,predict_grid_prototype,predict_cluster_prototype
      real(kind=wp),dimension(kohonen_map%number_variables1,kohonen_map%number_variables2) :: current_values
      real(kind=wp),dimension(kohonen_map%number_clusters) :: distance_units
   !
      number_variables=kohonen_map%number_variables1*kohonen_map%number_variables2;
   !
         write(*,*)
         write(*,*) ' TWO LEVEL SOM: Prediction starting...';
         write(*,*)
   !       write(*,*) '                Prediction for Grid Layer in progress'
         do ipattern = 1, size(input_data,1)
            ihit = 0;
            jhit = 0;
            khit = 0;
            dist_hit = 100000.0;
            call input_data(ipattern)%get(current_prototype);
            !call current_prototype%get_prototype(current_values);         
            do iz = 1, size(kohonen_map%grid,3);  
               do iy = 1, size(kohonen_map%grid,2);
                  do ix = 1,size(kohonen_map%grid,1);
                     dist = 0.0_wp;
                     dist=kohonen_map%grid(ix,iy,iz)%distance(current_prototype,kohonen_map%distance_function);
                     dist = dist/float(number_variables);
                     if (dist .lt. dist_hit) then
                        dist_hit = dist
                        ihit = ix;
                        jhit = iy;
                        khit = iz;
                     endif
                  enddo!ix
               enddo!iy
            enddo!iz
            !call kohonen_map%grid(ihit,jhit,khit)%get_prototype(current_values);
            predict_grid_prototype=kohonen_map%grid(ihit,jhit,khit);
            ihit=0;
            dist_hit=10.0e4;
            do ic=1,size(kohonen_map%cluster_layer)
               dist=kohonen_map%cluster_layer(ic)%distance(predict_grid_prototype,kohonen_map%distance_function);
               distance_units(ic)=dist;
               if(dist .lt. dist_hit) then
                  dist_hit=dist;
                  ihit=ic;
               endif            
            enddo
   !         write(*,*) ipattern,ihit
            predict_cluster_prototype=kohonen_map%cluster_layer(ihit);
            call predict_cluster_prototype%get_prototype(current_values);
            map_output(ipattern,1)=ihit
         enddo !ipattern
         write(*,*) 'TWO LEVEL SOM: Prediction finished';
   
   !
   end subroutine predict_2lsom